write.morpho: Generate a file with a morphological alignment for MCMCtree

View source: R/morpho.R

write.morphoR Documentation

Generate a file with a morphological alignment for MCMCtree

Description

Generate an alignment file with quantitative characters in MCMCtree format. The option "seqfile" in the control file used by MCMCtree should read the path to the file output by this function.

Usage

write.morpho(
  M,
  filename,
  c = 0,
  R = diag(1, dim(M)[2]),
  method,
  A = NULL,
  names,
  ages,
  scale = 1
)

Arguments

M

Matrix, s rows (specimens) and n morphological continuous characters (see details).

filename

Character, name for the output file.

c

Numeric, vector of variances within the species of a population (see details). If not provided, c = 0 (no population noise).

R

Matrix, correlation matrix. Requires method (see details). If not provided, R = I (no character correlation).

method

(Optional) character, either "eigen" or "chol", method used to decompose the inverse of the shrinkage correlation matrix. Requires R (see details).

A

(Optional) matrix, decomposed matrix. Requires R but not method (see details).

names

(Optional) list or character, species name included in the morphological alignment (see examples B and C).

ages

(Optional) list or numeric, ages of the species included in the morphological alignment (see example C).

scale

numeric, all characters are multiplied by scale before the morphological matrix is printed to file. Useful to re-scale the characters so they have the required variance (rate).

Details

The matrix M has s rows, one for each specimen, and n columns regarding the characters. If the data set contains landmarks, they 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 ...
... ... ... ... ... ... ... ...

See descriptions in data-raw/carnivores19x29.raw.R, data-raw/vulpes21x29.raw.R, data-raw/C.R, data-raw/V.R, to know how to obtain the morphological alignment used in write.morpho. Note that the explanation starts with the processing of raw data.

If the data set contains a set of n morphological continuous characters, e.g. from a simulated data set, the file should look like

specimens char.1 char.2 char.3 char.4 ...
Sp_1 0.143 -0.028 -0.044 0.129 ...
Sp_2 0.128 -0.024 -0.028 0.124 ...
... ... ... ... ... ...

Note that if a list with the specimens names is not passed to the parameter names, the name for each species will be "Species_1", "Species_2", and so on.

The object c can be of length 1, if all characters have the same variance, or a vector of length n with the variance of each of the characters. For the latter, you can take a look at object var.foxes, which has been generated following the steps explained in data-raw/var.foxes.R.

The object R has to be a symmetric and positive definite object of class matrix, i.e. class( R ) = "matrix". See R.sh for an example of its format. You can also read the description in data-raw/R.shrinkage.R for the details about how to generate this matrix.

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.

If a correlation matrix R is provided, write.morpho can use either the method = "chol" or method = "eigen" to get a matrix A such that \mathrm{R^{-1}}=\mathrm{A^{T}}\mathrm{A}. This matrix \mathrm{A^{T}} is later used to transform the morphological data to account for the correlation in this data set, so that the transformed characters in Z, \mathrm{Z}=\mathrm{M}\mathrm{A^{T}}, are independent. Alternatively, this matrix A can also be provided by the user. You can read the description to generate this matrix in data-raw/A.R to understand how to generate this matrix, which is also available as object A. If you decide to use a matrix A, it will be used to transform the data and no decomposition will be performed, thus saving computational time when large matrices are to be used.

Author(s)

Sandra Alvarez-Carretero and Mario dos Reis

See Also

matrix2array, array2matrix, sim.morpho

Examples

# A.1) Providing the morphological alignment (M = C) and
#      the name for the output file. This does not account for
#      correlation nor population noise

       write.morpho(  M = C, filename = "seqfile.aln" )

# A.2) Providing the morphological alignment (M = C), the population
#      noise (c = 0.25), and the name for the output file. Note that
#      c = 0.25 means that the population noise for all the characters
#      is c = 0.25, i.e. it will be considered as if
#      length( c ) = p characters, being all of them 0.25.

       write.morpho(  M = C, c = 0.25,
                      filename = "seqfile.aln" )

# A.3) Providing the morphological alignment (M = C), the population
#      noise (c = 0.25), the (estimate of the) correlation matrix (R),
#      the method to decompose R ("chol" in this example),
#      and the name for the output file. Note that the R matrix needs
#      to be invertible, otherwise the data will not be able to be
#      transformed accounting for correlation.

       write.morpho(  M = C, c = 0.25, R = R.sh,
                      method = "chol", filename = "seqfile.aln" )

# A.4) Providing the morphological alignment (M = C), a vector with
#      the population noise for each character (c = var.foxes),
#      the (shrinkage estimate of the) correlation matrix (R = R.sh),
#      the method to decompose R ("chol" in this example),
#      and the name for the output file. Note that the matrix passed
#      to argument R needs to be invertible, otherwise the data will
#      not be able to be transformed accounting for correlation.

       write.morpho(  M = C, c = var.foxes, R = R.sh,
                      method = "chol", filename = "seqfile.aln" )

# A.5) Providing the morphological alignment (C), a vector with
#      the population noise for each character (c = var.foxes),
#      the (shrinkage estimate of the) correlation matrix (R = R.sh),
#      the A matrix to transform the data, and the name for the
#      output file. Note that as the A matrix is provided, the matrix
#      passed to R will not be decomposed, hence the argument "method"
#      is not needed.

       write.morpho(  M = C, c = var.foxes, R = R.sh,
                      A = A, filename = "seqfile.aln" )

# B) Scenario A.5 but providing a list with the
#    names of the species

     names <- list( sp1  = "Ael_sp.", sp2  = "Can_dir", sp3  = "Epi_hay", sp4  = "Hes_sp.",
                    sp5  = "Par_jos", sp6  = "Tom_sp.", sp7  = "Enh_pah", sp8  = "Cuo_alp",
                    sp9  = "Spe_ven", sp10 = "Can_lup", sp11 = "Cer_tho", sp12 = "Oto_meg",
                    sp13 = "Urs_ame", sp14 = "Ail_ful", sp15 = "Nan_bin", sp16 = "Par_her",
                    sp17 = "Hia_won", sp18 = "Smi_fat", sp19 = "Vul_vul"
                  )

     write.morpho( M = C, c = var.foxes, R = R.sh,
                   A = A, filename = "seqfile.aln", names = names )

# C) Scenario A.5 but providing a vector of type character with the names of
#    the specimens and a list with their corresponding ages. Please
#    keep the same order in both lists, so the first specimen in the
#    list name corresponds to the first age in the age list, and so on.

     names <- c( "Ael_sp.", "Can_dir", "Epi_hay", "Hes_sp.",
                 "Par_jos", "Tom_sp.", "Enh_pah", "Cuo_alp",
                 "Spe_ven", "Can_lup", "Cer_tho", "Oto_meg",
                 "Urs_ame", "Ail_ful", "Nan_bin", "Par_her",
                 "Hia_won", "Smi_fat", "Vul_vul"
                  )

     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
                  )

     write.morpho( M = C, c = var.foxes, R = R.sh,
                   A = A, filename = "seqfile.aln",
                   names = names, ages = ages )


dosreislab/mcmc3r documentation built on July 21, 2024, 3:23 a.m.