BioExplorer.Community_MatrixType
mutable struct Community_Matrix

The Community_Matrix is a mutable composite object designed to represent data pertinent to taxonomic diversity (TD) analyses within the BioExplorer.jl package. It comprises four essential components:

Fields

  • sites::Vector{String}: A vector of strings representing the unique identifiers or names of ecological sites or locations where data on species are collected.
  • species::Vector{String}: A vector of strings denoting the unique identifiers or names of species observed or recorded within the ecological sites.
  • species_data::Array{Float32, 2}: An array of single-precision floating-point numbers (Float32) containing the actual data on species occurrences or abundances. This array is structured such that rows correspond to different ecological sites, while columns correspond to different species.
  • type::String: A string indicating the type of the data (either 'occurrences' or 'abundances') stored within the object.
BioExplorer.Trait_MatrixType
mutable struct Trait_Matrix

The Trait_Matrix is a mutable composite object developed for conducting functional diversity (FD) analyses within the BioExplorer.jl package. It encompasses four key components:

Fields

  • traits::Vector{String}: A vector of strings representing the unique identifiers or names of ecological traits under consideration for analysis.
  • species::Vector{String}: A vector of strings signifying the unique identifiers or names of species associated with the ecological traits.
  • species_data::Array{Any, 2}: An array of arbitrary data types (Any) containing the actual trait data corresponding to each species. This array is structured such that rows correspond to different traits, while columns represent different species.
  • type::Vector{String}: A vector of strings indicating the type of traits. Accepted values are "N", "C", or "O", respectively referring to Nominal or Binary traits, Continuous traits, and Ordinal traits.
BioExplorer.FD_dispersionMethod
FD_dispersion(trait_matrix::Trait_Matrix, weight)

Compute the functional dispersion of a trait matrix using Gower distance and PCA.

Arguments

  • trait_matrix::Trait_Matrix: A trait matrix containing trait values for different species.
  • weight: Vector of weight parameters for each traits in the trait matrix. Can be missing, in this case, each traits are consider equals.

Returns

  • dispersion: Functional dispersion value.

Details

The function computes the Gower distance matrix from the trait matrix and performs Principal Component Analysis (PCA) on it. It then projects the data onto the first two principal components and computes the convex hull of the projected points. Next, it calculates the centroid of the convex hull and computes the distance of each point on the hull from the centroid. The functional dispersion is then calculated as the mean distance of all points on the hull from the centroid.

BioExplorer.FD_obs_metricsMethod
FD_obs_metrics(trait_matrix::Trait_Matrix, weight)

Compute observation-based metrics of functional diversity: Originality, Uniqueness, and Contribution.

Arguments

  • trait_matrix::Trait_Matrix: A trait matrix containing trait values for different species.
  • weight: Vector of weight parameters for each traits in the trait matrix. Can be missing, in this case, each traits are consider equals.

Returns

  • Obs_metrics: DataFrame containing observation-based metrics for each species.

Details

The function calculates three observation-based metrics of functional diversity:

  1. Species Contribution: Contribution of each species to the overall functional richness.
  2. Species Originality: Average distance of each species to all other species in trait space.
  3. Species Uniqueness: Minimum distance of each species to any other species in trait space.

The contribution of each species is calculated by iteratively removing one species at a time, recalculating the functional richness, and computing the difference from the total richness. Originality is calculated as the average distance of each species to all other species, while uniqueness is calculated as the minimum distance of each species to any other species.

BioExplorer.FD_richMethod
FD_rich(trait_matrix::Trait_Matrix, weight, display_graph::Bool)

Calculate the functional richness of a trait matrix using Gower distance and PCA.

Arguments

  • trait_matrix::Trait_Matrix: A trait matrix containing trait values for different species.
  • weight: Vector of weight parameters for each traits in the trait matrix. Can be missing, in this case, each traits are consider equals.
  • display_graph::Bool: A boolean flag indicating whether to display the graph.

Returns

  • rich: Functional richness value.

Details

The function computes the Gower distance matrix from the trait matrix and performs Principal Component Analysis (PCA) on it. It then projects the data onto the first two principal components and computes the convex hull of the projected points. If display_graph is set to true, it generates a scatter plot of the projected points with the convex hull highlighted. Finally, it calculates the functional richness using Gauss' shoelace formula for the area of the convex hull.

BioExplorer.SACFunction
SAC(community_matrix::Community_Matrix, npermut::Int64 = 0)

Compute the species accumulation curve for a community matrix.

Arguments

  • community_matrix::Community_Matrix: The input community matrix.
  • npermut::Int64: The number of permutations for computing confidence intervals around the mean accumulation curve. Defaults to 0, indicating no permutation.

Returns

  • If npermut is 0, returns a vector representing the species accumulation curve.

  • If npermut is greater than 0, returns a DataFrame containing the mean accumulation curve and confidence intervals.

  • Makie Figure representing the species accumulation curve.

Details

The function computes the species accumulation curve for the input community matrix. If npermut is greater than 0, it performs permutations of the community matrix and computes multiple accumulation curves to estimate confidence intervals around the mean accumulation curve. The function generates a plot of the accumulation curve and displays it.

BioExplorer.beta_carvalhoMethod
beta_carvalho(community_matrix::Community_Matrix)

Compute beta diversity between two communities of a community matrix according to the Carvalho framework:

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species data for two communities. If more than communities are present, only the two first ones are considered.

Returns

  • An array containing three components of beta diversity: total dissimilarity, replacement dissimilarity, and richness dissimilarity.

Details

The function calculates beta diversity between two communities using the Carvalho framework, which decomposes beta diversity into three components: total dissimilarity, replacement dissimilarity, and richness dissimilarity. Total dissimilarity measures the overall difference in species composition between two communities. Replacement dissimilarity measures the extent to which species in one community are replaced by different species in the other community. Richness dissimilarity measures the difference in species richness between the two communities.

Reference

Carvalho, J. C., Cardoso, P., Borges, P. A. V., Schmera, D., & Podani, J. (2013). Measuring fractions of beta diversity and their relationships to nestedness: A theoretical and empirical comparison of novel approaches. Oikos, 122(6), 825-834. https://doi.org/10.1111/j.1600-0706.2012.20980.x

BioExplorer.beta_carvalho_matrixMethod
beta_carvalho_matrix(community_matrix::Community_Matrix)

Compute the beta diversity matrix between all pairs of communities in a community matrix according to the Carvalho framework: Carvalho, J. C., Cardoso, P., Borges, P. A. V., Schmera, D., & Podani, J. (2013). Measuring fractions of beta diversity and their relationships to nestedness: A theoretical and empirical comparison of novel approaches. Oikos, 122(6), 825-834. https://doi.org/10.1111/j.1600-0706.2012.20980.x

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species data for multiple communities.

Returns

  • An array containing three DataFrames representing the beta diversity matrices for total dissimilarity, replacement dissimilarity, and richness dissimilarity.

Details

The function calculates beta diversity between all pairs of communities using the Carvalho framework, which decomposes beta diversity into three components: total dissimilarity, replacement dissimilarity, and richness dissimilarity. For each component, a separate DataFrame is returned where each row and column represent a community, and the values represent the beta diversity between corresponding pairs of communities.

See also beta_carvalho

BioExplorer.biodiversity_estimateMethod

biodiversityestimate(communitymatrix::Community_Matrix)

Compute several estimations of the number of species in a community.

Arguments

  • community_matrix::Community_Matrix: The input community matrix.

Returns

  • A DataFrame containing the estimators and their corresponding values.

Details

This function computes several estimations of the number of species in a community based on incidence or abundance data. It calculates the Jackknife1, Jackknife2, Chao1, and Chao2 estimators of species richness.

See also sampling_coverage

Reference

Chao, A., Chiu, C.-H., & Jost, L. (2014). Unifying Species Diversity, Phylogenetic Diversity, Functional Diversity, and Related Similarity and Differentiation Measures Through Hill Numbers. Annual Review of Ecology, Evolution, and Systematics, 45(1), 297–324. https://doi.org/10.1146/annurev-ecolsys-120213-091540

BioExplorer.biodiversity_surfaceMethod
biodiversity_surface(community_matrix::Community_Matrix)

Create a biodiversity surface plot for a community matrix.

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species abundance data for multiple communities.

Returns

  • A Makie Figure representing the biodiversity surface plot.

Details

The function generates a 3D surface plot using Makie to visualize biodiversity across communities. It computes the Hill series for each community using the BioExplorer.hill function and plots the equivalent number of species against the order of Hill number. The x-axis represents different communities, the y-axis represents the order of Hill number, and the z-axis represents the equivalent number of species. Different colors on the surface indicate variations in biodiversity across communities.

See also hill

BioExplorer.fit_gambinMethod
fit_gambin(community_matrix::Community_Matrix, community_selected::String)

Fit a Gamma Binomial (GamBin) model to a species abundance distribution.

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species abundance data.
  • community_selected::String: The name of the community for which the model will be fitted.

Returns

  • The optimized parameter (α) of the Gamma Binomial model.

Details

The function fits a Gamma Binomial (GamBin) model to the species abundance distribution of the selected community. The model is fitted using maximum likelihood estimation to find the optimal value of the parameter (α) that maximizes the likelihood of observing the given abundance distribution under the GamBin model.

See also octave, octave_plot, fit_gambin_plot

Reference

Ugland, K. I., Lambshead, P. J. D., McGill, B., Gray, J. S., O’Dea, N., Ladle, R. J., & Whittaker, R. J. (2007). Modelling dimensionality in species abundance distributions: Description and evaluation of the Gambin model. Evolutionary Ecology Research, 9(2), 313–324.

Matthews, T. J., Borregaard, M. K., Ugland, K. I., Borges, P. A. V., Rigal, F., Cardoso, P., & Whittaker, R. J. (2014). The gambin model provides a superior fit to species abundance distributions with a single free parameter: Evidence, implementation and interpretation. Ecography, 37(10), 1002–1011. https://doi.org/10.1111/ecog.00861

BioExplorer.fit_gambin_plotMethod
fit_gambin_plot(community_matrix::Community_Matrix, community_selected::String)

Plot the species abundance distribution and the fitted Gamma Binomial (GamBin) model.

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species abundance data.
  • community_selected::String: The name of the community for which the model will be fitted and plotted.

Returns

  • A Makie Figure representing the plot of species abundance distribution and the fitted GamBin model.

Details

The function plots the species abundance distribution of the selected community along with the fitted Gamma Binomial (GamBin) model. It first fits the GamBin model using the fit_gambin function to obtain the optimized parameter (α), and then plots the observed abundance distribution as a bar plot and the fitted model as a line plot on the same axis. The resulting Makie Figure represents the plot of species abundance distribution and the fitted GamBin model.

See also octave, octave_plot, fit_gambin

Reference

Ugland, K. I., Lambshead, P. J. D., McGill, B., Gray, J. S., O’Dea, N., Ladle, R. J., & Whittaker, R. J. (2007). Modelling dimensionality in species abundance distributions: Description and evaluation of the Gambin model. Evolutionary Ecology Research, 9(2), 313–324.

Matthews, T. J., Borregaard, M. K., Ugland, K. I., Borges, P. A. V., Rigal, F., Cardoso, P., & Whittaker, R. J. (2014). The gambin model provides a superior fit to species abundance distributions with a single free parameter: Evidence, implementation and interpretation. Ecography, 37(10), 1002–1011. https://doi.org/10.1111/ecog.00861

BioExplorer.generate_abundance_communitiesFunction
generate_abundance_communities(nb_species::Int64, nb_sites::Int64, poisson_lambda::Number = Float64(50))

Generate a fake community matrix with abundance data using Poisson distribution.

Arguments

  • nb_species::Int64: The number of species in each community.
  • nb_sites::Int64: The number of sites (communities).
  • poisson_lambda::Number: The lambda parameter for the Poisson distribution used to generate abundance data. Defaults to 50.

Returns

  • A community matrix containing randomly generated abundance data.

Details

The function generates a fake community matrix with abundance data by randomly sampling species abundances from a Poisson distribution. Each row of the matrix represents a site (community), and each column represents a species. The abundance data is rounded to the nearest integer. The default value of the lambda parameter for the Poisson distribution is set to 50.

See also generate_incidence_communities

BioExplorer.generate_incidence_communitiesFunction
generate_incidence_communities(nb_species::Int64, nb_sites::Int64, probability::Number = Float64(50))

Generate a fake community matrix with incidence data using a Binomial distribution.

Arguments

  • nb_species::Int64: The number of species in each community.
  • nb_sites::Int64: The number of sites (communities).
  • probability::Number: The probability parameter for the Binomial distribution used to generate incidence data. Defaults to 50.

Returns

  • A community matrix containing randomly generated incidence data.

Details

The function generates a fake community matrix with incidence data by randomly sampling species presence/absence using a Binomial distribution. Each row of the matrix represents a site (community), and each column represents a species. The incidence data is rounded to the nearest integer (0 or 1). The default value of the probability parameter for the Binomial distribution is set to p = 0.5.

See also generate_abundance_communities

BioExplorer.hillMethod
hill(community_matrix::Community_Matrix)

Compute the Hill series for all communities within a community matrix.

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species abundance data for multiple communities.

Returns

  • A DataFrame representing the Hill series for each community, including Hill numbers H0, H1, H2, and H3.

Details

The function calculates the Hill series, which represents the equivalent number of species at different orders of Hill number, for each community within the given community matrix. Hill numbers measure biodiversity by scaling species abundance to a common unit, providing a way to compare biodiversity across different communities. The output DataFrame contains the Hill numbers H0, H1, H2, and H3 for each community.

Reference

Hill, M. O. (1973). Diversity and Evenness: A Unifying Notation and Its Consequences. Ecology, 54(2), 427–432. https://doi.org/10.2307/1934352

BioExplorer.jaccard_dissimMethod
jaccard_dissim(community_matrix::Community_Matrix)

Compute the Jaccard dissimilarity between two communities.

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species data for two communities. If more than communities are present, only the two first ones are considered.

Returns

  • jaccard_dissim_value::Float64: The Jaccard dissimilarity value between the two communities.

Details

The Jaccard dissimilarity measures dissimilarity between two sets by comparing their intersection to their union. For two communities represented by binary species presence-absence data, it is calculated as the number of species found in only one of the communities divided by the total number of species found in either community.

BioExplorer.jaccard_dissim_matrixMethod
jaccard_dissim_matrix(community_matrix::Community_Matrix)

Compute the Jaccard dissimilarity matrix between all communities in a community matrix.

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species data for multiple communities.

Returns

  • beta_matrix::DataFrame: A DataFrame representing the Jaccard dissimilarity matrix between all pairs of communities.

Details

The function calculates the Jaccard dissimilarity between all pairs of communities in the given community matrix. It iterates over all possible combinations of communities and computes the Jaccard dissimilarity using the jaccard_dissim function. The result is stored in a DataFrame where each row and column represent a community, and the values represent the Jaccard dissimilarity between corresponding pairs of communities.

See also jaccard_dissim

BioExplorer.mat_com_convertMethod
mat_com_convert(community_matrix::Community_Matrix)

Convert a community matrix with abundance data to a community matrix with incidence data.

Arguments

  • community_matrix::Community_Matrix: The input community matrix with abundance data.

Returns

  • The corresponding community matrix with incidence data.

Details

The function converts a community matrix with abundance data to a community matrix with incidence data by replacing abundance values greater than 0 with 1 (indicating presence) and leaving abundance values of 0 unchanged (indicating absence). The resulting matrix has the same structure as the input matrix but with incidence data.

BioExplorer.matrix_GowdisMethod
matrix_Gowdis(trait_matrix::Trait_Matrix, weight)

Compute the Gower distance between every pair of species within a trait matrix.

Arguments

  • trait_matrix::Trait_Matrix: A trait matrix containing species trait data.
  • weight: Vector of weight parameters for each traits in the trait matrix. Can be missing, in this case, each traits are consider equals.

Returns

  • A matrix representing the Gower distance between every pair of species.

Details

The function computes the Gower distance between every pair of species based on their trait data within the given trait matrix. It iterates over all combinations of species and computes the pairwise Gower distance using the pairwise_Gowdis function. The resulting matrix represents the Gower distance between every pair of species.

See also pairwise_Gowdis

BioExplorer.octaveMethod
octave(community_matrix::Community_Matrix, community_selected::String)

Compute the octave distribution of species in a selected community.

Arguments

  • community_matrix::Community_Matrix: The input community matrix.
  • community_selected::String: The name of the community for which the octave distribution is to be computed.

Returns

A tuple (community_name, community_ranked) where:

  • community_name is the name of the selected community.
  • community_octave is a DataFrame containing the species abundance and their corresponding octave class.

Details

This function computes the octave distribution of species abundance within a selected community from the input community matrix. The octave distribution categorizes species abundance into octave classes, with each class representing a doubling in abundance.

See also octave_plot,

BioExplorer.octave_plotMethod
octave_plot(community_matrix::Community_Matrix, community_selected::String)

Display the Octave plot of a selected community.

Arguments

  • community_matrix::Community_Matrix: The input community matrix.
  • community_selected::String: The name of the community for which the Octave plot is to be displayed.

Returns

  • Makie Figure representing the Octave plot.

Details

This function computes the octave distribution of species abundance within a selected community using the octave function and displays it as a bar plot. Each bar in the plot represents an octave class, with the x-axis indicating the octave class and the y-axis indicating the number of species falling within each octave class.

See also octave,

BioExplorer.pairwise_GowdisMethod
pairwise_Gowdis(trait_matrix::Trait_Matrix, weight, species1::String, species2::String)

Compute the Gower distance between a pair of species within a trait matrix.

Arguments

  • trait_matrix::Trait_Matrix: A trait matrix containing species trait data.
  • weight: Vector of weight parameters for each traits in the trait matrix. Can be missing, in this case, each traits are consider equals.
  • species1::String: The name of the first species.
  • species2::String: The name of the second species.

Returns

  • The Gower distance between the two species.

Details

The function computes the Gower distance between a pair of species based on their trait data within the given trait matrix. Gower distance is a measure of dissimilarity that accounts for different types of variables (categorical, numerical, and ordinal) and their respective weights. The dissimilarity is calculated for each trait variable and then combined using the provided weights to obtain the final Gower distance.

See also matrix_Gowdis

BioExplorer.pielou_evennessMethod
pielou_evenness(community_matrix::Community_Matrix)

Compute an evenness diversity metric according to the Pielou's framework.

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species abundance data.

Returns

  • A DataFrame containing the Pielou's evenness (J) metric for each community.

Details

The function calculates the Pielou's evenness (J) metric for each community within the given community matrix. Pielou's evenness measures the equitability of species abundance distribution within a community. It is calculated as the Shannon-Wiener entropy divided by the maximum possible entropy, which is the logarithm of the species richness.

Reference

Pielou, E. C. (1969). An Introduction to Mathematical Ecology. Wiley-Interscience.

BioExplorer.rankMethod
rank(community_matrix::Community_Matrix, community_selected::String)

Rank the species according to their abundance in a selected community.

Arguments

  • community_matrix::Community_Matrix: The input community matrix.
  • community_selected::String: The name of the community for which species ranking is required.

Returns

A tuple (community_name, community_ranked) where:

  • community_name is the name of the selected community.
  • community_ranked is a DataFrame containing species ranked by abundance in the selected community.

Details

This function ranks the species according to their abundance in the specified community of the input community matrix. It sorts the species by abundance in descending order and assigns ranks. If species have the same abundance, they are assigned the same rank (ex-aequo ranking). The function removes species with zero abundance.

See also whittacker_plot

BioExplorer.sampling_coverageMethod
sampling_coverage(community_matrix::Community_Matrix)

Compute the percentage of sampling coverage according to the Chao framework.

Arguments

  • community_matrix::Community_Matrix: A community matrix containing species abundance data.

Returns

  • The percentage of sampling coverage.

Details

The function calculates the percentage of sampling coverage, which estimates the proportion of total species in the community that have been sampled. It uses the Chao framework, which considers the abundance of rare species to estimate the coverage. The output represents the percentage of species that have been sampled based on the observed data.

Reference

Chao, A., Kubota, Y., Zelený, D., Chiu, C., Li, C., Kusumoto, B., Yasuhara, M., Thorn, S., Wei, C., Costello, M. J., & Colwell, R. K. (2020). Quantifying sample completeness and comparing diversities among assemblages. Ecological Research, 35(2), 292–314. https://doi.org/10.1111/1440-1703.12102

BioExplorer.whittacker_plotMethod
whittacker_plot(community_matrix::Community_Matrix, community_selected::String)

Display the Whittaker plot of a selected community.

Arguments

  • community_matrix::Community_Matrix: The input community matrix.
  • community_selected::String: The name of the community for which the Whittaker plot is to be displayed.

Returns

  • Makie Figure representing the Whittaker plot.

Details

This function generates a Whittaker plot for the selected community from the input community matrix. The plot shows the log10 abundance of species against their rank. The species are ranked by abundance, and the abundance is plotted on the y-axis in log10 scale. The rank of each species is plotted on the x-axis.

See also [rank](@re, octave_plot