proc2MCMCtree: Procrustes alignment output in MCMCtree format for Bayesian...

View source: R/morpho.R

proc2MCMCtreeR Documentation

Procrustes alignment output in MCMCtree format for Bayesian inference

Description

Perform a Procrustes alignment with the procSym and align2procSym and generate a morphological alignment in MCMCtree format.

Usage

proc2MCMCtree(
  data,
  popdata,
  sp.data,
  sp.popdata,
  filename,
  coords,
  method = c("eigen", "chol"),
  ages,
  ...
)

Arguments

data

Matrix or array, object with the data set with one specimen per species. See details.

popdata

Matrix or array, containing the specimens for the population of species with a larger variation (more than one specimen for this species). See details.

sp.data

Numeric, position of the specimen in the data object also present in the popdata object (see details).

sp.popdata

Numeric, position of the specimen in the popdata object also present in the data object (see details).

filename

Character, name for the output file.

coords

Numeric, 2 or 3 for 2D or 3D landmarks, respectively.

method

(Optional) character, either "eigen" or "col", method used to decompose the inverse of the shrinkage correlation matrix. See details.

ages

(Optional) list or vector, ages of the species included in the morpholical alignment.

...

(Optional) Further arguments passed to procSym (see details)

Details

If the objects data and popdata are of class "matrix", then they have s or ps rows, respectively, and n columns as of the amount of characters. These data sets are supposed to contain landmarks, which can be given in 2D or 3D. For instance, if the landmarks are 3D, the first 3 columns will be the coordinates x, y, and z for the first landmark, the next 3 columns for the second landmark; and so on.

specimens lmk1.x lmk1.y lmk1.z lmk2.x lmk2.y lmk2.z ...
Sp_1 0.143 -0.028 -0.044 0.129 0.028 -0.043 ...
Sp_2 0.128 -0.024 -0.028 0.124 0.027 -0.025 ...
... ... ... ... ... ... ... ...

You can load the object C.mat.unal to see an example of this format and also take a look at the data-raw/C.R file for details about how to generate it. Otherwise, if the objects data and popdata are of class "array", they have format k x q x s or k x q x ps, respectively, where k is the number of landmarks, q the number of coordinates, s the number of specimens for object data, and ps the number of specimens in a sampled population of one species for object popdata. You can load the object C.arr.unal to see an example of this array format and also take a look at the data-raw/C.R file for details about how to generate it.

The names of the specimens will be taken from either the row names of the data matrices, if data and popdata objects are of class "matrix", or from the names of the third dimension of the array, if these objects are of class "array". If no names are found, the specimens will be labelled as "Specimen_1", "Specimen_2", and so on.

Note that you are providing a population matrix with the landmarks collected from more than on specimen belonging to the same species. Therefore, it is assumed that population noise is accounted for. Furthermore, the landmarks of one of these specimens are also present in the morphological alignment. First, a Procrusts alignment is performed with the data set with one specimen per species (data), and then all the specimens in popdata except for the specimen common in data, which position is specified in the argument sp.popdata, are aligned to the mean shape of the PA generated with data. Take a look at data.raw/V.R for more details.

The logarithm of the determinant of the correlation matrix is going to be printed in the output file to later be used by MCMCtree during the likelihood calculation.

The function procSym can perform a Procrustes superimposition alignment given a dataset of landmarks in array format and allows the user to pass different options to the arguments defined. The current function runs with the default parameters in procSym but allows the user to pass three arguments to proc2MCMCtree through ...: scale, reflect, and pairedLM. If you are thinking of using these arguments, read the documentation in procSym, otherwise do not pass further arguments to proc2MCMCtree and let it run procSym in default mode.

Value

$dataPS

Object of class symproc output by the function Morpho::procSym. The user can access the rest of the items in this object by checking <your_object>$dataPS$<available_Morpho_objects>. It contains the array with the morphological alignment generated by this function, which corresponds to the specimens in the object passed to data. This can be accessed by loading <your_object>$dataPS$rotated.

$popdataPS

Array with the morphological alignment with the object passed to popdata, i.e. the alignment with the population sample

$M

Matrix with the morphological alignment with the specimens of all species. Note that this alignment is not corrected for character correlation nor population noise. The corrected alignment is only printed in the output alignment file.

$Rsh

Estimated shrinkage correlation matrix.

$c

Estimated population variance.

Author(s)

Sandra Alvarez-Carretero and Mario dos Reis

See Also

write.morpho, matrix2array, array2matrix

Examples

## Not run: 
# A. Use the unaligned, but processed, carnivoran data (data = C.mat.unal) and
#    vulpes data (popdata = V.mat.unal) to obtain a morphological alignment.
#    The fox specimen that is common in V.mat.unal and C.mat.unal is the
#    one in the first row of V.mat.unal (sp.popdata = 1) and the one in position
#    13 in C.mat.unal (sp.data = 13). The method to
#    decompose the estimate shrinkage correlation matrix (internally calculated
#    in this function) is the Cholesky decomposition (method = c("chol")).
#    We do not add ages.
     right    <- c( 11, 22, 13, 19, 15, 20, 24, 5, 7, 2, 9, 26, 29 )
     left     <- c( 10, 21, 12, 17, 14, 18, 23, 4, 6, 1, 8, 25, 28 )
     pairedLM <- cbind( left, right )
     obj.aln <- proc2MCMCtree( data = C.mat.unal, popdata = V.mat.unal, sp.data = 13,
                               sp.popdata = 1, filename = "./seqfile.aln", coords = 3,
                               method = c("chol"), pairedLM = pairedLM )

# B. Use the unaligned, but processed, carnivoran data (data = C.mat.unal) and
#    vulpes data (popdata = V.mat.unal) to obtain a morphological alignment.
#    The fox specimen that is common in V.mat.unal and C.mat.unal is the
#    one in the first row of V.mat.unal (sp.popdata = 1) and the one in position
#    13 in C.mat.unal (sp.data = 13). The method to
#    decompose the estimate shrinkage correlation matrix (internally calculated
#    in this function) is the Cholesky decomposition (method = c("chol")).
#    We add ages.
     ages <- list( sp1  = 13.135, sp2  =  0.0285, sp3  = 11.95,  sp4  = 35.55,
                   sp5  = 25.615, sp6  = 14.785,  sp7  = 28.55,  sp8  =  0,
                   sp9  =  0,     sp10 =  0,      sp11 =  0,     sp12 =  0,
                   sp13 =  0,     sp14 =  0,      sp15 =  0,     sp16 =  0,
                   sp17 =  6.65,  sp18 =  0.0285, sp19 =  0 )
     right    <- c( 11, 22, 13, 19, 15, 20, 24, 5, 7, 2, 9, 26, 29 )
     left     <- c( 10, 21, 12, 17, 14, 18, 23, 4, 6, 1, 8, 25, 28 )
     pairedLM <- cbind( left, right )
     obj.aln <- proc2MCMCtree( data = C.mat.unal, popdata = V.mat.unal, sp.data = 13,
                               sp.popdata = 1, filename = "./seqfile.aln", coords = 3,
                               method = c("chol"), pairedLM = pairedLM, ages = ages )

## End(Not run)


dosreislab/mcmc3r documentation built on March 29, 2024, 6:45 p.m.